Contracting cubeΒΆ
In this demo we simulate a unit cube that is fixed at \(x = 0\) and free at \(x = 1\). We use a transversally isotropic material with fiber oriented in the \(x\)-direction.
import dolfin
try:
from dolfin_adjoint import (
Constant,
DirichletBC,
Expression,
UnitCubeMesh,
interpolate,
Mesh,
)
except ImportError:
from dolfin import (
Constant,
DirichletBC,
interpolate,
Expression,
UnitCubeMesh,
Mesh,
)
import pulse
from fenics_plotly import plot
pulse.iterate.logger.setLevel(10)
# Create mesh
N = 6
mesh = UnitCubeMesh(N, N, N)
# Create subdomains
class Free(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] > (1.0 - dolfin.DOLFIN_EPS) and on_boundary
class Fixed(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] < dolfin.DOLFIN_EPS and on_boundary
# Create a facet fuction in order to mark the subdomains
ffun = dolfin.MeshFunction("size_t", mesh, 2)
ffun.set_all(0)
# Mark the first subdomain with value 1
fixed = Fixed()
fixed_marker = 1
fixed.mark(ffun, fixed_marker)
# Mark the second subdomain with value 2
free = Free()
free_marker = 2
free.mark(ffun, free_marker)
# Create a cell function (but we are not using it)
cfun = dolfin.MeshFunction("size_t", mesh, 3)
cfun.set_all(0)
# Collect the functions containing the markers
marker_functions = pulse.MarkerFunctions(ffun=ffun, cfun=cfun)
# Create mictrotructure
V_f = dolfin.VectorFunctionSpace(mesh, "CG", 1)
# Fibers
f0 = interpolate(Expression(("1.0", "0.0", "0.0"), degree=1), V_f)
# Sheets
s0 = interpolate(Expression(("0.0", "1.0", "0.0"), degree=1), V_f)
# Fiber-sheet normal
n0 = interpolate(Expression(("0.0", "0.0", "1.0"), degree=1), V_f)
# Collect the mictrotructure
microstructure = pulse.Microstructure(f0=f0, s0=s0, n0=n0)
# Create the geometry
geometry = pulse.Geometry(
mesh=mesh,
marker_functions=marker_functions,
microstructure=microstructure,
)
# Use the default material parameters
material_parameters = pulse.HolzapfelOgden.default_parameters()
# Select model for active contraction
active_model = pulse.ActiveModels.active_strain
# active_model = "active_stress"
# Set the activation
activation = Constant(0.0)
# Create material
material = pulse.HolzapfelOgden(
active_model=active_model,
parameters=material_parameters,
activation=activation,
)
# Make Dirichlet boundary conditions
def dirichlet_bc(W):
V = W if W.sub(0).num_sub_spaces() == 0 else W.sub(0)
return DirichletBC(V, Constant((0.0, 0.0, 0.0)), fixed)
# Make Neumann boundary conditions
neumann_bc = pulse.NeumannBC(traction=Constant(0.0), marker=free_marker)
# Collect Boundary Conditions
bcs = pulse.BoundaryConditions(dirichlet=(dirichlet_bc,), neumann=(neumann_bc,))
# Create problem
problem = pulse.MechanicsProblem(geometry, material, bcs)
# Solve problem
pulse.iterate.iterate(problem, activation, 0.1)
2021-11-27 09:59:41,793 [995] DEBUG pulse.iterate: Control: [0.0]
2021-11-27 09:59:41,794 [995] DEBUG pulse.iterate: Target: [0.1]
2021-11-27 09:59:41,796 [995] DEBUG pulse.iterate: Intial number of steps: 5 with step size 0.02
2021-11-27 09:59:41,796 [995] INFO pulse.iterate: Iterating to target control (f_21)...
2021-11-27 09:59:41,797 [995] INFO pulse.iterate: Current control: f_21 = 0.000
2021-11-27 09:59:41,797 [995] INFO pulse.iterate: Target: 0.100
2021-11-27 09:59:41,798 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 09:59:41,799 [995] DEBUG pulse.iterate: Maximum difference: 1.000e-01
2021-11-27 09:59:41,799 [995] DEBUG pulse.iterate: Try new control
2021-11-27 09:59:41,800 [995] DEBUG pulse.iterate: Current control: 0.020
2021-11-27 10:00:42,695 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:00:42,697 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:42,698 [995] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2021-11-27 10:00:42,701 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:42,702 [995] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2021-11-27 10:00:42,703 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:42,704 [995] DEBUG pulse.iterate: Current control: 0.040
2021-11-27 10:00:47,036 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:00:47,037 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:47,037 [995] DEBUG pulse.iterate: Maximum difference: 6.000e-02
2021-11-27 10:00:47,039 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:47,039 [995] DEBUG pulse.iterate: Maximum difference: 6.000e-02
2021-11-27 10:00:47,040 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:47,042 [995] DEBUG pulse.iterate: Current control: 0.060
2021-11-27 10:00:49,135 [995] DEBUG pulse.iterate: Solver did not converge...
2021-11-27 10:00:49,136 [995] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-27 10:00:49,137 [995] DEBUG pulse.iterate: Reduce control step
2021-11-27 10:00:49,138 [995] DEBUG pulse.iterate: Assign old state
2021-11-27 10:00:49,146 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:49,147 [995] DEBUG pulse.iterate: Maximum difference: 6.000e-02
2021-11-27 10:00:49,149 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:49,150 [995] DEBUG pulse.iterate: Current control: 0.050
2021-11-27 10:00:51,678 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:00:51,679 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:51,682 [995] DEBUG pulse.iterate: Maximum difference: 5.000e-02
2021-11-27 10:00:51,683 [995] DEBUG pulse.iterate: Adapt step size. New step size: 0.015
2021-11-27 10:00:51,685 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:51,685 [995] DEBUG pulse.iterate: Maximum difference: 5.000e-02
2021-11-27 10:00:51,687 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:51,688 [995] DEBUG pulse.iterate: Current control: 0.065
2021-11-27 10:00:53,877 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:00:53,878 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:53,879 [995] DEBUG pulse.iterate: Maximum difference: 3.500e-02
2021-11-27 10:00:53,880 [995] DEBUG pulse.iterate: Adapt step size. New step size: 0.022
2021-11-27 10:00:53,881 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:53,882 [995] DEBUG pulse.iterate: Maximum difference: 3.500e-02
2021-11-27 10:00:53,883 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:53,883 [995] DEBUG pulse.iterate: Current control: 0.087
2021-11-27 10:00:55,509 [995] DEBUG pulse.iterate: Solver did not converge...
2021-11-27 10:00:55,510 [995] DEBUG pulse.iterate:
NOT CONVERGING
2021-11-27 10:00:55,511 [995] DEBUG pulse.iterate: Reduce control step
2021-11-27 10:00:55,512 [995] DEBUG pulse.iterate: Assign old state
2021-11-27 10:00:55,520 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:55,520 [995] DEBUG pulse.iterate: Maximum difference: 3.500e-02
2021-11-27 10:00:55,521 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:55,522 [995] DEBUG pulse.iterate: Current control: 0.076
2021-11-27 10:00:57,491 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:00:57,492 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:57,493 [995] DEBUG pulse.iterate: Maximum difference: 2.375e-02
2021-11-27 10:00:57,494 [995] DEBUG pulse.iterate: Adapt step size. New step size: 0.017
2021-11-27 10:00:57,495 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:57,496 [995] DEBUG pulse.iterate: Maximum difference: 2.375e-02
2021-11-27 10:00:57,498 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:57,498 [995] DEBUG pulse.iterate: Current control: 0.093
2021-11-27 10:00:59,657 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:00:59,658 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:59,660 [995] DEBUG pulse.iterate: Maximum difference: 6.875e-03
2021-11-27 10:00:59,661 [995] DEBUG pulse.iterate: Adapt step size. New step size: 0.025
2021-11-27 10:00:59,662 [995] DEBUG pulse.iterate: Check target reached: NO
2021-11-27 10:00:59,663 [995] DEBUG pulse.iterate: Maximum difference: 6.875e-03
2021-11-27 10:00:59,664 [995] DEBUG pulse.iterate: Change step size for final iteration
2021-11-27 10:00:59,665 [995] DEBUG pulse.iterate: Try new control
2021-11-27 10:00:59,666 [995] DEBUG pulse.iterate: Current control: 0.100
2021-11-27 10:01:01,266 [995] DEBUG pulse.iterate:
SUCCESFULL STEP:
2021-11-27 10:01:01,267 [995] DEBUG pulse.iterate: Check target reached: YES!
2021-11-27 10:01:01,269 [995] DEBUG pulse.iterate: Check target reached: YES!
([Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 47),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 152),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 227),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 321),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 364),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 442),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 485),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 520)],
[Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 45),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 150),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 225),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 319),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 362),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 440),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 483),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 518)])
# Get displacement and hydrostatic pressure
u, p = problem.state.split(deepcopy=True)
V = dolfin.VectorFunctionSpace(mesh, "CG", 1)
u_int = interpolate(u, V)
new_mesh = Mesh(mesh)
dolfin.ALE.move(new_mesh, u_int)
fig = plot(mesh, show=False)
fig.add_plot(plot(new_mesh, color="red", show=False))
fig.show()